Field Sources

The superscreen.sources module provides factory functions for a few convenient superscreen.parameter.Parameter classes used to define applied magnetic fields.

ConstantField

superscreen.sources.ConstantField(value=0)[source]

Returns a Parameter that computes a constant as a function of x, y, z.

Parameters:

value (float) – The constant value of the field.

Return type:

Parameter

Returns:

A Parameter that returns value at all x, y, z.

MonopoleField

superscreen.sources.MonopoleField(r0=(0, 0, 0), nPhi0=1, vector=False)[source]

Returns a Parameter that computes the z-component of the field from a monopole (monopole) located at position (x0, y0, z0) containing a total of nPhi0 flux quanta.

\[\mu_0H_z(\vec{r}-\vec{r}_0) = \frac{n\Phi_0}{2\pi} \frac{(\vec{r}-\vec{r}_0)\cdot\hat{z}}{|(\vec{r}-\vec{r}_0)|^3}\]
Parameters:
  • r0 (Tuple[float, float, float]) – Coordinates of the monopole position.

  • nPhi0 (Union[int, float]) – Number of flux quanta contained in the monopole.

  • vector (bool) – If True, return the vector magnetic field. Otherwise, return only the \(z\)-component.

Return type:

Parameter

Returns:

A Parameter that returns the field in units of Phi_0 / (length_units)**2.

PearlVortexField

superscreen.sources.PearlVortexField(*, r0=(0, 0, 0), Lambda=0, nPhi0=1, xs, ys)[source]

Returns a Parameter that computes the z-component of the field from a Pearl vortex located at position (x0, y0, z0) in a film with effective penetration depth Lambda (Pearl length 2 * Lambda) containing a total of nPhi0 flux quanta.

The field from a Pearl vortex located at is computed using a Fourier transform method. For a uniform thin film lying in the \(x-y\) plane with effective penetration depth \(\Lambda\) (Pearl length \(2\Lambda\)), the Fourier transform of the \(z\)-component of the field from a vortex containing \(n\) flux quanta located at the origin, \(x=y=z=0\), is given by:

\[\mathcal{F}\{\mu_0H_z\}(k_x, k_y, z) = \frac{n\Phi_0e^{-kz}}{1 + 2\Lambda k},\]

where \(k=\sqrt{k_x^2 + k_y^2}\) and the quantity is in units of Phi_0 / (length_units)**2, where length_units are the units of xs, ys, etc. The field is calculated by inverse Fourier-transforming the above expression for an \(x-y\) plane defined by parameters xs and ys, then interpolating the field to the desired coordinates. Note that the Fourier method may not be accurate if xs and ys are not sampled finely enough.

See also

References: [Pearl-APL-1964], [Tafuri-PRL-2004].

Parameters:
  • r0 (Tuple[float, float, float]) – Coordinates of the Pearl vortex position.

  • Lambda (float) – The effective penetration depth of the film in which the vortex lies. Lambda is equal to half the Pearl length.

  • nPhi0 (Union[int, float]) – Number of flux quanta contained in the monopole.

  • xs (ndarray) – Vectors of x and y coordinates defining the the domain in which the field will be computed using a Fourier transform as described above.

  • ys (ndarray) – Vectors of x and y coordinates defining the the domain in which the field will be computed using a Fourier transform as described above.

Return type:

Parameter

Returns:

A Parameter that returns the out-of-plane field in units of Phi_0 / (length_units)**2.

DipoleField

superscreen.sources.DipoleField(*, dipole_positions, dipole_moments, component=None, length_units='um', moment_units='mu_B')[source]

Returns a Parameter that computes a given component of the field from a distribution of dipoles with given moments (in units of the Bohr magneton) located at the given positions.

Given dipole positions \(\vec{r}_{0, i}\) and moments \(\vec{m}_i\), the magnetic field is:

\[\mu_0\vec{H}(\vec{r}) = \sum_i\frac{\mu_0}{4\pi} \frac{ 3\hat{r}_i(\hat{r}_i\cdot\vec{m}_i) - \vec{m}_i }{ |\vec{r}_i|^3 },\]

where \(\vec{r}_i=(x, y, z) - \vec{r}_{0, i}\).

Parameters:
  • dipole_positions (Union[ndarray, Tuple[float, float, float]]) – Coordinates (x0_i, y0_i, z0_i) of the position of each dipole i. Shape (3, ) or (1, 3) for a single dipole, or shape (m, 3) for m dipoles.

  • dipole_moments (Union[ndarray, Tuple[float, float, float]]) – Dipole moments (mx_i, my_i, mz_i) in units of the Bohr magneton. If dipole_moments has shape (3, ) or (1, 3), then all dipoles are assigned the same moment. Otherwise, dipole_moments must have shape (m, 3), i.e. the moment is specified for each dipole.

  • component (Optional[str]) – The component of the field to calculate: “x”, “y”, “z”, or None. If None, then the vector field (shape (m, 3)) is returned.

  • length_units (str) – The units for the positions coordinates x, y, z, and dipole_positions.

  • moment_units (str) – The units for dipole_moments, for example the Bohr magneton “mu_B” or SI base units “A * m ** 2”.

Return type:

Parameter

Returns:

A Parameter that computes a given component of the field \(\mu_0\vec{H}(x, y, z)\) in Tesla for a given distribution of dipoles.

superscreen.sources.dipole_field(eval_coords, r0=(0, 0, 0), moment=(0, 0, 0))[source]

Returns the 3D field from a single dipole with the given moment (in units of amps * meters ** 2) located at the position r0, evaluated at coordinates eval_coords = [x, y, z].

Given \(\vec{r}=(x, y, z) - \vec{r}_0\), the magnetic field is:

\[\vec{B}(\vec{r}) = \frac{\mu_0}{4\pi} \frac{ 3\hat{r}(\hat{r}\cdot\vec{m}) - \vec{m} }{ |\vec{r}|^3 }\]
Parameters:
  • eval_coords (ndarray) – (x, y, z) coordinates (in meters) at which to evaluate the field. Either a sequence of length 3 (for a single position) or an array of shape (n, 3) (for n positions.).

  • r0 (Union[ndarray, Tuple[float, float, float]]) – Coordinates (x0, y0, z0) (in meters) of the dipole position, shape (3,) or (1, 3).

  • moment (Union[ndarray, Tuple[float, float, float]]) – Dipole moment (mx, my, mz) in units of amps * meters ** 2, shape (3,) or (1, 3).

Returns:

An array with shape (3, ) if x, y, z are scalars, or shape (n, 3) if x, y, z are vectors with shape (n, ).

Return type:

Magnetic field (Bx, By, Bz) in Tesla evaluated at (x, y, z)

superscreen.sources.dipole_distribution(x, y, z, *, dipole_positions, dipole_moments, component=None, length_units='um', moment_units='mu_B')[source]

Returns the 3D field \(\vec{B}=\mu_0\vec{H}\), or one of its components, from a distribution of dipoles with given moments (in units of the Bohr magneton) located at the given positions, evaluated at coordinates (x, y, z).

Parameters:
  • x (Union[float, ndarray]) – x-coordinate(s) (in meters) at which to evaluate the field. Either a scalar or vector with shape (n, ).

  • y (Union[float, ndarray]) – y-coordinate(s) (in meters) at which to evaluate the field. Either a scalar or vector with shape (n, ).

  • z (Union[float, ndarray]) – z-coordinate(s) (in meters) at which to evaluate the field. Either a scalar or vector with shape (n, ).

  • dipole_positions (ndarray) – Coordinates (x0_i, y0_i, z0_i) of the position of each dipole i, shape (m, 3) (in meters) .

  • dipole_moments (Union[ndarray, Tuple[float, float, float]]) – Dipole moments (mx_i, my_i, mz_i) in units of amps * meters ** 2. If dipole_moments has shape (3, ) or (1, 3), then all dipoles are assigned the same moment. Otherwise, dipole_moments must have shape (m, 3), i.e. the moment is specified for each dipole.

  • component (Optional[str]) – The component of the magnetic field to return: “x”, “y”, “z”, or None. If None, the vector magnetic field (shape (n, 3)) is returned.

  • length_units (str) – The units for the positions coordinates x, y, z, and dipole_positions.

  • moment_units (str) – The units for dipole_moments, for example the Bohr magneton “mu_B” or SI base units “A * m ** 2”.

Return type:

ndarray

Returns:

Magnetic field (Bx, By, Bz) (or one of its components) in Tesla evaluated at (x, y, z): An array with shape (3, ) if x, y, z are scalars, or shape (n, 3) if x, y, z are vectors with shape (n, ).

SheetCurrentField

superscreen.sources.SheetCurrentField(*, sheet_positions, current_densities, z0, length_units='um', current_units='uA')[source]

Returns a Parameter that computes the z-component of the field from a 2D sheet of current parameterized by the given positions and current densities.

The \(z\)-component of the field from a 2D sheet of current \(S\) lying in the plane \(z=z_0\) with spatially varying current density \(\vec{J}=(J_x, J_y)\) is given by:

\[\mu_0H_z(\vec{r})=\frac{\mu_0}{4\pi}\int_S \frac{J_x(\vec{r}')(\vec{r}-\vec{r}')\cdot\hat{y} - J_y(\vec{r}')(\vec{r}-\vec{r}')\cdot\hat{x}} {|\vec{r}-\vec{r}'|^3}\,\mathrm{d}^2r',\]

where \(\vec{r}=(x, y, z)\) and \(\vec{r}'=(x', y', z_0)\).

Parameters:
  • sheet_positions (ndarray) – Coordinates (x0, y0) (in meters) of the current sheet, shape (m, 2).

  • current_densities (ndarray) – 2D current density (Jx, Jy) in units of amps / meter, shape (m, 2).

  • z0 (float) – Vertical (z) position of the current sheet.

  • length_units (str) – The units for all coordinates.

  • current_units (str) – The units for current values. The current_densities are assumed to be in units of current_units / length_units.

Return type:

Parameter

Returns:

A Parameter that computes \(\mu_0\vec{H}_z(x, y, z)\) in Tesla for a given sheet current.

superscreen.sources.biot_savart_2d(x, y, z, *, positions, current_densities, z0=0, areas=None, length_units='um', current_units='uA', vector=True)[source]

Returns the magnetic field (in tesla) from a sheet of current located at vertical positon z0 (in units of length_units). The current is parameterized by a set of current_densities (in units of current_units / length_units) and x-y positions (in units of length_units), and the field is evaluated at coordinates (x, y, z).

\[\begin{split}\mu_0H_x(\vec{r}) &= \frac{\mu_0}{4\pi}\int_S \frac{J_y(\vec{r}')(\vec{r}-\vec{r}')\cdot\hat{z}} {|\vec{r}-\vec{r}'|^3}\,\mathrm{d}^2r'\\ \mu_0H_y(\vec{r}) &= \frac{\mu_0}{4\pi}\int_S -\frac{J_x(\vec{r}')(\vec{r}-\vec{r}')\cdot\hat{z}} {|\vec{r}-\vec{r}'|^3}\,\mathrm{d}^2r'\\ \mu_0H_z(\vec{r}) &= \frac{\mu_0}{4\pi}\int_S \frac{J_x(\vec{r}')(\vec{r}-\vec{r}')\cdot\hat{y} - J_y(\vec{r}')(\vec{r}-\vec{r}')\cdot\hat{x}} {|\vec{r}-\vec{r}'|^3}\,\mathrm{d}^2r'\end{split}\]

where \(\vec{r}=(x, y, z)\) and \(\vec{r}'=(x', y', z_0)\).

Parameters:
  • x (Union[float, ndarray]) – x-coordinate(s) at which to evaluate the field. Either a scalar or vector with shape (n, ).

  • y (Union[float, ndarray]) – y-coordinate(s) at which to evaluate the field. Either a scalar or vector with shape (n, ).

  • z (Union[float, ndarray]) – z-coordinate(s) at which to evaluate the field. Either a scalar or vector with shape (n, ).

  • positions (ndarray) – Coordinates (x0, y0) of the current sheet, shape (m, 2).

  • current_densities (ndarray) – 2D current density (Jx, Jy), shape``(m, 2)``.

  • z0 (float) – Vertical (z) position of the current sheet.

  • areas (Optional[ndarray]) – Vertex areas for positions in units of length_units**2. If None, the positions are triangulated to calculate vertex areas.

  • length_units (str) – The units for all coordinates.

  • current_units (str) – The units for current values. The current_densities are assumed to be in units of current_units / length_units.

  • vector (bool) – Return the full vector magnetic field (shape (n, 3)) rather than just the z-component (shape (n, )).

Return type:

ndarray

Returns:

Magnetic field in tesla evaluated at (x, y, z). If vector is True, returns the vector magnetic field \(\mu_0\vec{H}\) (shape (n, 3)). Otherwise, returns the the \(z\)-component, \(\mu_0H_z\) (shape (n,)).