simsopt.field package

Submodules

simsopt.field.biotsavart module

class simsopt.field.biotsavart.BiotSavart(coils)

Bases: BiotSavart, MagneticField

Computes the MagneticField induced by a list of closed curves \(\Gamma_k\) with electric currents \(I_k\). The field is given by

\[B(\mathbf{x}) = \frac{\mu_0}{4\pi} \sum_{k=1}^{n_\mathrm{coils}} I_k \int_0^1 \frac{(\Gamma_k(\phi)-\mathbf{x})\times \Gamma_k'(\phi)}{\|\Gamma_k(\phi)-\mathbf{x}\|^3} d\phi\]

where \(\mu_0=4\pi 10^{-7}\) is the magnetic constant.

Parameters:

coils – A list of simsopt.field.coil.Coil objects.

A_and_dA_vjp(v, vgrad)

Same as simsopt.geo.biotsavart.BiotSavart.A_vjp but returns the vector Jacobian product for \(A\) and \(\nabla A\), i.e. it returns

\[\{ \sum_{i=1}^{n} \mathbf{v}_i \cdot \partial_{\mathbf{c}_k} \mathbf{A}_i \}_k, \{ \sum_{i=1}^{n} {\mathbf{v}_\mathrm{grad}}_i \cdot \partial_{\mathbf{c}_k} \nabla \mathbf{A}_i \}_k.\]
A_vjp(v)

Assume the field was evaluated at points \(\mathbf{x}_i, i\in \{1, \ldots, n\}\) and denote the value of the field at those points by \(\{\mathbf{A}_i\}_{i=1}^n\). These values depend on the shape of the coils, i.e. on the dofs \(\mathbf{c}_k\) of each coil. This function returns the vector Jacobian product of this dependency, i.e.

\[\{ \sum_{i=1}^{n} \mathbf{v}_i \cdot \partial_{\mathbf{c}_k} \mathbf{A}_i \}_k.\]
B_and_dB_vjp(v, vgrad)

Same as simsopt.geo.biotsavart.BiotSavart.B_vjp but returns the vector Jacobian product for \(B\) and \(\nabla B\), i.e. it returns

\[\{ \sum_{i=1}^{n} \mathbf{v}_i \cdot \partial_{\mathbf{c}_k} \mathbf{B}_i \}_k, \{ \sum_{i=1}^{n} {\mathbf{v}_\mathrm{grad}}_i \cdot \partial_{\mathbf{c}_k} \nabla \mathbf{B}_i \}_k.\]
B_vjp(v)

Assume the field was evaluated at points \(\mathbf{x}_i, i\in \{1, \ldots, n\}\) and denote the value of the field at those points by \(\{\mathbf{B}_i\}_{i=1}^n\). These values depend on the shape of the coils, i.e. on the dofs \(\mathbf{c}_k\) of each coil. This function returns the vector Jacobian product of this dependency, i.e.

\[\{ \sum_{i=1}^{n} \mathbf{v}_i \cdot \partial_{\mathbf{c}_k} \mathbf{B}_i \}_k.\]
as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

d2A_by_dXdcoilcurrents(compute_derivatives=1)
d2B_by_dXdcoilcurrents(compute_derivatives=1)
d3A_by_dXdXdcoilcurrents(compute_derivatives=2)
d3B_by_dXdXdcoilcurrents(compute_derivatives=2)
dA_by_dcoilcurrents(compute_derivatives=0)
dB_by_dcoilcurrents(compute_derivatives=0)
classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

simsopt.field.boozermagneticfield module

class simsopt.field.boozermagneticfield.BoozerAnalytic(etabar, B0, N, G0, psi0, iota0, Bbar=1.0, I0=0.0, G1=0.0, I1=0.0, K1=0.0)

Bases: BoozerMagneticField

Computes a BoozerMagneticField based on a first-order expansion in distance from the magnetic axis (Landreman & Sengupta, Journal of Plasma Physics 2018). Here the magnetic field strength is expressed as,

\[B(s,\theta,\zeta) = B_0 \left(1 + \overline{\eta} \sqrt{2s\psi_0/\overline{B}}\cos(\theta - N \zeta)\right),\]

the covariant components are,

\[ \begin{align}\begin{aligned}G(s) = G_0 + \sqrt{2s\psi_0/\overline{B}} G_1\\I(s) = I_0 + \sqrt{2s\psi_0/\overline{B}} I_1\\K(s,\theta,\zeta) = \sqrt{2s\psi_0/\overline{B}} K_1 \sin(\theta - N \zeta),\end{aligned}\end{align} \]

and the rotational transform is,

\[\iota(s) = \iota_0.\]

While formally \(I_0 = I_1 = G_1 = K_1 = 0\), these terms have been included in order to test the guiding center equations at finite beta.

Parameters:
  • etabar – magnitude of first order correction to magnetic field strength

  • B0 – magnetic field strength on the axis

  • N – helicity of symmetry (integer)

  • G0 – lowest order toroidal covariant component

  • psi0 – (toroidal flux)/ (2*pi) on the boundary

  • iota0 – lowest order rotational transform

  • Bbar – normalizing magnetic field strength (defaults to 1)

  • I0 – lowest order poloidal covariant component (defaults to 0)

  • G1 – first order correction to toroidal covariant component (defaults to 0)

  • I1 – first order correction to poloidal covariant component (defaults to 0)

  • K1 – first order correction to radial covariant component (defaults to 0)

set_B0(B0)
set_Bbar(Bbar)
set_G0(G0)
set_G1(G1)
set_I0(I0)
set_I1(I1)
set_K1(K1)
set_N(N)
set_etabar(etabar)
set_iota0(iota0)
set_psi0(psi0)
class simsopt.field.boozermagneticfield.BoozerMagneticField(psi0)

Bases: BoozerMagneticField

Generic class that represents a magnetic field in Boozer coordinates \((s,\theta,\zeta)\). Here \(s = \psi/\psi_0\) is the normalized toroidal flux where \(2\pi\psi_0\) is the toroidal flux at the boundary. The magnetic field in the covariant form is,

\[\textbf B(s,\theta,\zeta) = G(s) \nabla \zeta + I(s) \nabla \theta + K(s,\theta,\zeta) \nabla \psi,\]

and the contravariant form is,

\[\textbf B(s,\theta,\zeta) = \frac{1}{\sqrt{g}} \left(\frac{\partial \mathbf r}{\partial \zeta} + \iota(s)\frac{\partial \mathbf r}{\partial \theta}\right),\]

where,

\[\sqrt{g}(s,\theta,\zeta) = \frac{G(s) + \iota(s)I(s)}{B^2}.\]

Here \(\iota(s) = \psi_P'(\psi)\) where \(2\pi\psi_P\) is the poloidal flux and \(2\pi\psi\) is the toroidal flux. Each subclass of BoozerMagneticField implements functions to compute \(B\), \(G\), \(I\), \(\iota\), \(\psi_P\), and their derivatives. The cylindrical coordinates \(R(s,\theta,\zeta)\) and \(Z(s,\theta,\zeta)\) in addition to \(K(s,\theta,\zeta)\) and \(\nu\) where \(\zeta = \phi + \nu(s,\theta,\zeta)\) and \(\phi\) is the cylindrical azimuthal angle are also implemented by BoozerRadialInterpolant and InterpolatedBoozerField. The usage is similar to the MagneticField class.

The usage of BoozerMagneticField` is as follows:

booz = BoozerAnalytic(etabar,B0,N,G0,psi0,iota0) # An instance of BoozerMagneticField
points = ... # points is a (n, 3) numpy array defining :math:`(s,\theta,\zeta)`
booz.set_points(points)
modB = bfield.modB() # returns the magnetic field strength at `points`

BoozerMagneticField has a cache to avoid repeated calculations. To clear this cache manually, call the clear_cached_properties() function. The cache is automatically cleared when set_points() is called or one of the dependencies changes.

clear_cached_properties()

Clear the cache.

recompute_bell(parent=None)
class simsopt.field.boozermagneticfield.BoozerRadialInterpolant(equil, order, mpol=32, ntor=32, N=None, enforce_vacuum=False, rescale=False, ns_delete=0, no_K=False)

Bases: BoozerMagneticField

Given a Vmec instance, performs a Boozer coordinate transformation using BOOZXFORM. The magnetic field can be computed at any point in Boozer coordinates using radial spline interpolation (scipy.interpolate.InterpolatedUnivariateSpline) and an inverse Fourier transform in the two angles. Throughout stellarator symmetry is assumed.

Parameters:
  • equil – instance of simsopt.mhd.vmec.Vmec or simsopt.mhd.boozer.Boozer. If it is an instance of simsopt.mhd.boozer.Boozer, the compute_surfs needs to include all of the grid points in the half-radius grid of the corresponding Vmec equilibrium.

  • order – (int) order for radial interpolation. Must satisfy 1 <= order <= 5.

  • mpol – (int) number of poloidal mode numbers for BOOZXFORM (defaults to 32)

  • ntor – (int) number of toroidal mode numbers for BOOZXFORM (defaults to 32)

  • N – Helicity of quasisymmetry to enforce. If specified, then the non-symmetric Fourier harmonics of \(B\) and \(K\) are filtered out. Otherwise, all harmonics are kept. (defaults to None)

  • enforce_vacuum – If True, a vacuum field is assumed, \(G\) is set to its mean value, \(I = 0\), and \(K = 0\).

  • rescale

    If True, use the interpolation method in the DELTA5D code. Here, a few of the first radial grid points or (bmnc, rmnc, zmns, numns, kmns) are deleted (determined by ns_delete). The Fourier harmonics are then rescaled as:

    bmnc(s)/s^(1/2) for m = 1

    bmnc(s)/s for m even and >= 2

    bmnc(s)/s^(3/2) for m odd and >=3

    before performing interpolation and spline differentiation to obtain dbmncds. If False, interpolation of the unscaled Fourier harmonics and its finite-difference derivative wrt s is performed instead (defaults to False)

  • ns_delete – (see rescale) (defaults to 0)

compute_K()
init_splines()
class simsopt.field.boozermagneticfield.InterpolatedBoozerField(field, degree, srange, thetarange, zetarange, extrapolate=True, nfp=1, stellsym=True)

Bases: InterpolatedBoozerField, BoozerMagneticField

This field takes an existing BoozerMagneticField and interpolates it on a regular grid in \(s,\theta,\zeta\). This resulting interpolant can then be evaluated very quickly. This is modeled after InterpolatedField.

__init__(field, degree, srange, thetarange, zetarange, extrapolate=True, nfp=1, stellsym=True)
Parameters:
  • field – the underlying simsopt.field.boozermagneticfield.BoozerMagneticField to be interpolated.

  • degree – the degree of the piecewise polynomial interpolant.

  • srange – a 3-tuple of the form (smin, smax, ns). This mean that the interval [smin, smax] is split into ns many subintervals.

  • thetarange – a 3-tuple of the form (thetamin, thetamax, ntheta). thetamin must be >= 0, and thetamax must be <=2*pi.

  • zetarange – a 3-tuple of the form (zetamin, zetamax, nzeta). zetamin must be >= 0, and thetamax must be <=2*pi.

  • extrapolate – whether to extrapolate the field when evaluate outside the integration domain or to throw an error.

  • nfp – Whether to exploit rotational symmetry. In this case any toroidal angle is always mapped into the interval \([0, 2\pi/\mathrm{nfp})\), hence it makes sense to use zetamin=0 and zetamax=2*np.pi/nfp.

  • stellsym – Whether to exploit stellarator symmetry. In this case theta is always mapped to the interval \([0, \pi]\), hence it makes sense to use thetamin=0 and thetamax=np.pi.

simsopt.field.coil module

class simsopt.field.coil.Coil(curve, current)

Bases: Coil, Optimizable

A Coil combines a Curve and a Current and is used as input for a BiotSavart field.

plot(**kwargs)

Plot the coil’s curve. This method is just shorthand for calling the plot() function on the underlying Curve. All arguments are passed to simsopt.geo.curve.Curve.plot()

vjp(v_gamma, v_gammadash, v_current)
class simsopt.field.coil.Current(current, dofs=None, **kwargs)

Bases: Current, CurrentBase

An optimizable object that wraps around a single scalar degree of freedom. It represents the electric current in a coil, or in a set of coils that are constrained to use the same current.

property current
vjp(v_current)
simsopt.field.coil.apply_symmetries_to_currents(base_currents, nfp, stellsym)

Take a list of n Current`s and return ``n * nfp * (1+int(stellsym))` Current objects obtained by copying (for nfp rotations) and sign-flipping (optionally for stellarator symmetry).

simsopt.field.coil.apply_symmetries_to_curves(base_curves, nfp, stellsym)

Take a list of n simsopt.geo.curve.Curve`s and return ``n * nfp * (1+int(stellsym))` simsopt.geo.curve.Curve objects obtained by applying rotations and flipping corresponding to nfp fold rotational symmetry and optionally stellarator symmetry.

simsopt.field.coil.coils_to_focus(filename, curves, currents, nfp=1, stellsym=False, Ifree=False, Lfree=False)

Export a list of Curve objects together with currents in FOCUS format, so they can be used by FOCUS. The format is introduced at https://princetonuniversity.github.io/FOCUS/rdcoils.pdf This routine only works with curves of type CurveXYZFourier, not other curve types.

Parameters:
  • filename – Name of the file to write.

  • curves – A python list of CurveXYZFourier objects.

  • currents – Coil current of each curve.

  • nfp – The number of field periodicity. Defaults to 1.

  • stellsym – Whether or not following stellarator symmetry. Defaults to False.

  • Ifree – Flag specifying whether the coil current is free. Defaults to False.

  • Lfree – Flag specifying whether the coil geometry is free. Defaults to False.

simsopt.field.coil.coils_to_makegrid(filename, curves, currents, groups=None, nfp=1, stellsym=False)

Export a list of Curve objects together with currents in MAKEGRID input format, so they can be used by MAKEGRID and FOCUS. The format is introduced at https://princetonuniversity.github.io/STELLOPT/MAKEGRID Note that this function does not generate files with MAKEGRID’s output format.

Parameters:
  • filename – Name of the file to write.

  • curves – A python list of Curve objects.

  • currents – Coil current of each curve.

  • groups – Coil current group. Coils in the same group will be assembled together. Defaults to None.

  • nfp – The number of field periodicity. Defaults to 1.

  • stellsym – Whether or not following stellarator symmetry. Defaults to False.

simsopt.field.coil.coils_via_symmetries(curves, currents, nfp, stellsym)

Take a list of n curves and return n * nfp * (1+int(stellsym)) Coil objects obtained by applying rotations and flipping corresponding to nfp fold rotational symmetry and optionally stellarator symmetry.

simsopt.field.coil.load_coils_from_makegrid_file(filename, order, ppp=20)

This function loads a file in MAKEGRID input format containing the Cartesian coordinates and the currents for several coils and returns an array with the corresponding coils. The format is described at https://princetonuniversity.github.io/STELLOPT/MAKEGRID

Parameters:
  • filename – file to load.

  • order – maximum mode number in the Fourier expansion.

  • ppp – points-per-period: number of quadrature points per period.

Returns:

A list of Coil objects with the Fourier coefficients and currents given by the file.

simsopt.field.magneticfield module

class simsopt.field.magneticfield.MagneticField(**kwargs)

Bases: MagneticField, Optimizable

Generic class that represents any magnetic field from which each magnetic field class inherits. The usage of MagneticField is as follows:

bfield = BiotSavart(coils) # An instance of a MagneticField
points = ... # points is a (n, 3) numpy array
bfield.set_points(points)
B = bfield.B() # returns the Magnetic field at `points`
dA = bfield.dA_by_dX() # returns the gradient of the potential of the field at `points`

MagneticField has a cache to avoid repeated calculations. To clear this cache manually, call the clear_cached_properties() function. The cache is automatically cleared when set_points is called or one of the dependencies changes.

__add__(other)

Add two magnetic fields.

__mul__(other)

Multiply a field with a scalar.

__rmul__(other)

Multiply a field with a scalar.

clear_cached_properties()

Clear the cache.

recompute_bell(parent=None)

Function to be called whenever new DOFs input is given or if the parent Optimizable’s data changed, so the output from the current Optimizable object is invalid.

This method gets called by various DOF setters. If only the local DOFs of an object are being set, the recompute_bell method is called in that object and also in the descendent objects that have a dependency on the object, whose local DOFs are being changed. If gloabl DOFs of an object are being set, the recompute_bell method is called in the object, ancestors of the object, as well as the descendents of the object.

Need to be implemented by classes that provide a dof_setter for external handling of DOFs.

set_points(self: simsoptpp.MagneticField, arg0: numpy.ndarray[numpy.float64]) simsoptpp.MagneticField

Shorthand for set_points_cart.

set_points_cart(self: simsoptpp.MagneticField, arg0: numpy.ndarray[numpy.float64]) simsoptpp.MagneticField

Set the points where to evaluate the magnetic fields, in cartesian coordinates.

set_points_cyl(self: simsoptpp.MagneticField, arg0: numpy.ndarray[numpy.float64]) simsoptpp.MagneticField

Set the points where to evaluate the magnetic fields, in cylindrical coordinates (the order is \((r, \phi, z)\)).

to_vtk(filename, nr=10, nphi=10, nz=10, rmin=1.0, rmax=2.0, zmin=-0.5, zmax=0.5)

Export the field evaluated on a regular grid for visualisation with e.g. Paraview.

class simsopt.field.magneticfield.MagneticFieldMultiply(scalar, Bfield)

Bases: MagneticField

Class used to multiply a magnetic field by a scalar. It takes as input a MagneticField class and a scalar and multiplies B, A and their derivatives by that value.

as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

class simsopt.field.magneticfield.MagneticFieldSum(Bfields)

Bases: MagneticField

Class used to sum two or more magnetic field together. It can either be called directly with a list of magnetic fields given as input and outputing another magnetic field with B, A and its derivatives added together or it can be called by summing magnetic fields classes as Bfield1 + Bfield1

B_vjp(v)
as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

simsopt.field.magneticfieldclasses module

class simsopt.field.magneticfieldclasses.CircularCoil(r0=0.1, center=[0, 0, 0], I=159154.94309189534, normal=[0, 0])

Bases: MagneticField

Magnetic field created by a single circular coil evaluated using analytical functions, including complete elliptic integrals of the first and second kind. As inputs, it takes the radius of the coil (r0), its center, current (I) and its normal vector [either spherical angle components (normal=[theta,phi]) or (x,y,z) components of a vector (normal=[x,y,z])]). The (theta,phi) angles are related to the (x,y,z) components of the normal vector via theta = np.arctan2(normal[1], normal[0]) and phi = np.arctan2(np.sqrt(normal[0]**2+normal[1]**2), normal[2]). Sign convention: CircularCoil with a positive current produces a magnetic field vector in the same direction as the normal when evaluated at the center of the coil.a

Parameters:
  • r0 – radius of the coil

  • center – point at the coil center

  • I – current of the coil in Ampere’s

  • normal – if list with two values treats it as spherical angles theta and phi of the normal vector to the plane of the coil centered at the coil center, if list with three values treats it a vector

property I
as_dict(serial_objs_dict)

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

class simsopt.field.magneticfieldclasses.DipoleField(dipole_grid, dipole_vectors, stellsym=True, nfp=1, coordinate_flag='cartesian', m_maxima=None, R0=1)

Bases: MagneticField

Computes the MagneticField induced by N dipoles. The field is given by

\[B(\mathbf{x}) = \frac{\mu_0}{4\pi} \sum_{i=1}^{N} \left(\frac{3\mathbf{r}_i\cdot \mathbf{m}_i}{|\mathbf{r}_i|^5}\mathbf{r}_i - \frac{\mathbf{m}_i}{|\mathbf{r}_i|^3}\right)\]

where \(\mu_0=4\pi\times 10^{-7}\;N/A^2\) is the permeability of free space and \(\mathbf{r_i} = \mathbf{x} - \mathbf{x}^{dipole}_i\) is the vector between the field evaluation point and the dipole \(i\) position.

Parameters:
  • dipole_grid – 2D numpy array, shape (ndipoles, 3). A set of points corresponding to the locations of magnetic dipoles.

  • dipole_vectors – 2D numpy array, shape (ndipoles, 3). The dipole vectors of each of the dipoles in the grid.

  • stellsym – bool (default True). Whether or not the dipole grid is stellarator symmetric.

  • nfp – int (default 1). The field-period symmetry of the dipole-grid.

  • coordinate_flag – string (default “cartesian”). The global coordinate system that should be considered grid-aligned in the calculation. The options are “cartesian” (rectangular bricks), “cylindrical” (cylindrical bricks), and “toroidal” (uniform grid in simple toroidal coordinates). Note that this ASSUMES that the global coordinate system for the dipole locations is one of these three choices, so be careful if your permanent magnets are shaped/arranged differently!

  • m_maxima – 1D numpy array, shape (ndipoles,). The maximum dipole strengths of each magnet in the grid. If not specified, defaults to using the largest dipole strength of the magnets in dipole_grid, and using this value for all the dipoles. Needed for plotting normalized dipole magnitudes in the vtk functionality.

  • R0 – double. The value of the major radius of the stellarator needed only for simple toroidal coordinates.

_dipole_fields_from_symmetries(dipole_grid, dipole_vectors, stellsym=True, nfp=1, coordinate_flag='cartesian', m_maxima=None, R0=1)

Takes the dipoles and grid initialized in a PermanentMagnetOptimizer (for a half-period surface) and generates the full dipole manifold so that the call to B() (the magnetic field from the dipoles) correctly returns contributions from all the dipoles from symmetries.

_toVTK(vtkname)

Write dipole data into a VTK file (acknowledgements to Caoxiang’s CoilPy code).

Parameters:

vtkname (str) – VTK filename, will be appended with .vts or .vtu.

class simsopt.field.magneticfieldclasses.Dommaschk(mn=[[0, 0]], coeffs=[[0, 0]])

Bases: MagneticField

Vacuum magnetic field created by an explicit representation of the magnetic field scalar potential as proposed by W. Dommaschk (1986), Computer Physics Communications 40, 203-218. As inputs, it takes the arrays for the harmonics m, n and its corresponding coefficients.

Parameters:
  • m – first harmonic array

  • n – second harmonic array

  • coeffs – coefficient for Vml for each of the ith index of the harmonics m and n

as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

property mn
class simsopt.field.magneticfieldclasses.InterpolatedField(field, degree, rrange, phirange, zrange, extrapolate=True, nfp=1, stellsym=False, skip=None)

Bases: InterpolatedField, MagneticField

This field takes an existing field and interpolates it on a regular grid in \(r,\phi,z\). This resulting interpolant can then be evaluated very quickly.

__init__(field, degree, rrange, phirange, zrange, extrapolate=True, nfp=1, stellsym=False, skip=None)
Parameters:
  • field – the underlying simsopt.field.magneticfield.MagneticField to be interpolated.

  • degree – the degree of the piecewise polynomial interpolant.

  • rrange – a 3-tuple of the form (rmin, rmax, nr). This mean that the interval \([rmin, rmax]\) is split into nr many subintervals.

  • phirange – a 3-tuple of the form (phimin, phimax, nphi).

  • zrange – a 3-tuple of the form (zmin, zmax, nz).

  • extrapolate – whether to extrapolate the field when evaluate outside the integration domain or to throw an error.

  • nfp – Whether to exploit rotational symmetry. In this case any angle is always mapped into the interval \([0, 2\pi/\mathrm{nfp})\), hence it makes sense to use phimin=0 and phimax=2*np.pi/nfp.

  • stellsym – Whether to exploit stellarator symmetry. In this case z is always mapped to be positive, hence it makes sense to use zmin=0.

  • skip

    a function that takes in a point (in cylindrical (r,phi,z) coordinates) and returns whether to skip that location when building the interpolant or not. The signature should be

    def skip(r: double, phi: double, z: double) -> bool:
        ...
    

    See also here https://github.com/hiddenSymmetries/simsopt/pull/227 for a graphical illustration.

to_vtk(filename)

Export the field evaluated on a regular grid for visualisation with e.g. Paraview.

class simsopt.field.magneticfieldclasses.MirrorModel(B0=6.51292, gamma=0.124904, Z_m=0.98)

Bases: MagneticField

Model magnetic field employed in https://arxiv.org/abs/2305.06372 to study the magnetic mirror experiment WHAM. The magnetic field is given by \(\vec{B}=B_R \vec{e}_R + B_Z \vec{e}_Z\), where \(\vec{e}_R\) and \(\vec{e}_Z\) are the cylindrical radial and axial unit vectors, respectively, and \(B_R\) and \(B_Z\) are given by

\[B_R = -\frac{1}{R} \frac{\partial\psi}{\partial Z}, \; B_Z = \frac{1}{R}. \frac{\partial\psi}{\partial R}\]

In this model, the magnetic flux function \(\psi\) is written as a double Lorentzian function

\[\psi = \frac{R^2 \mathcal{B}}{2 \pi \gamma}\left(\left[1+\left(\frac{Z-Z_m}{\gamma}\right)^2\right]^{-1}+\left[1+\left(\frac{Z+Z_m}{\gamma}\right)^2\right]^{-1}\right).\]

Note that this field is neither a vacuum field nor a solution of MHD force balance. The input parameters are B0, gamma and Z_m with the standard values the ones used in https://arxiv.org/abs/2305.06372, that is, B0 = 6.51292, gamma = 0.124904, and Z_m = 0.98.

Parameters:
  • B0 – parameter \(\mathcal{B}\) of the flux surface function

  • gamma – parameter \(\gamma\) of the flux surface function

  • Z_m – parameter \(Z_m\) of the flux surface function

as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

class simsopt.field.magneticfieldclasses.PoloidalField(R0, B0, q)

Bases: MagneticField

Magnetic field purely in the poloidal direction, that is, in the theta direction of a poloidal-toroidal coordinate system. Its modulus is given by B = B0 * r / (R0 * q) so that, together with the toroidal field, it creates a safety factor equals to q

Parameters:
  • B0 – modulus of the magnetic field at R0

  • R0 – major radius of the magnetic axis

  • q – safety factor/pitch angle of the magnetic field lines

as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

class simsopt.field.magneticfieldclasses.Reiman(iota0=0.15, iota1=0.38, k=[6], epsilonk=[0.01], m0=1)

Bases: MagneticField

Magnetic field model in section 5 of Reiman and Greenside, Computer Physics Communications 43 (1986) 157—167. This field allows for an analytical expression of the magnetic island width that can be used for island optimization. However, the field is not completely physical as it does not have nested flux surfaces.

Parameters:
  • iota0 – unperturbed rotational transform

  • iota1 – unperturbed global magnetic shear

  • k – integer array specifying the Fourier modes used

  • epsilonk – coefficient of the Fourier modes

  • m0 – toroidal symmetry parameter (normally m0=1)

as_dict(serial_objs_dict)

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

class simsopt.field.magneticfieldclasses.ScalarPotentialRZMagneticField(phi_str)

Bases: MagneticField

Vacuum magnetic field as a solution of B = grad(Phi) where Phi is the magnetic field scalar potential. It takes Phi as an input string, which should contain an expression involving the standard cylindrical coordinates (R, phi, Z) Example: ScalarPotentialRZMagneticField(“2*phi”) yields a magnetic field B = grad(2*phi) = (0,2/R,0). In order for the analytical derivatives to be performed by sympy, a term 1e-30*Phi*R*Z is added to every entry. Note: this function needs sympy.

Parameters:

phi_str – string containing vacuum scalar potential expression as a function of R, Z and phi

as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

class simsopt.field.magneticfieldclasses.ToroidalField(R0, B0)

Bases: MagneticField

Magnetic field purely in the toroidal direction, that is, in the phi direction with (R,phi,Z) the standard cylindrical coordinates. Its modulus is given by B = B0*R0/R where R0 is the first input and B0 the second input to the function.

Parameters:
  • B0 – modulus of the magnetic field at R0

  • R0 – radius of normalization

as_dict(serial_objs_dict) dict

A JSON serializable dict representation of an object.

classmethod from_dict(d, serial_objs_dict, recon_objs)
Parameters:

d – Dict representation.

Returns:

GSONable class.

simsopt.field.normal_field module

class simsopt.field.normal_field.NormalField(nfp=1, stellsym=True, mpol=1, ntor=0, vns=None, vnc=None)

Bases: Optimizable

NormalField represents the magnetic field normal to a toroidal surface, for example the computational boundary of SPEC free-boundary.

Parameters:
  • nfp – The number of field period

  • stellsym – Whether (=True) or not (=False) stellarator symmetry is enforced.

  • mpol – Poloidal Fourier resolution

  • ntor – Toroidal Fourier resolution

  • vns

    Odd fourier modes of \(\mathbf{B}\cdot\mathbf{\hat{n}}\). 2D array of size (mpol+1)x(2ntor+1). Set to None to fill with zeros

    vns( mm, self.ntor+nn ) is the mode (mm,nn)

  • vnc

    Even fourier modes of \(\mathbf{B}\cdot\mathbf{\hat{n}}\). 2D array of size (mpol+1)x(2ntor+1). Ignored if stellsym if True. Set to None to fill with zeros

    vnc( mm, self.ntor+nn ) is the mode (mm,nn)

_make_names()

Form a list of names of the vns, vnc

change_resolution(mpol, ntor)

Change the values of mpol and ntor. Any new Fourier amplitudes will have a magnitude of zero. Any previous nonzero Fourier amplitudes that are not within the new range will be discarded.

check_mn(m, n)
fixed_range(mmin, mmax, nmin, nmax, fixed=True)

Set the ‘fixed’ property for a range of m and n values.

All modes with m in the interval [mmin, mmax] and n in the interval [nmin, nmax] will have their fixed property set to the value of the fixed parameter. Note that mmax and nmax are included (unlike the upper bound in python’s range(min, max).)

In case of non stellarator symmetric field, both vns and vnc will be fixed (or unfixed)

classmethod from_spec(filename)

Initialize using the harmonics in SPEC input file

get_index_in_dofs(m, n, mpol=None, ntor=None, even=False)

Returns position of mode (m,n) in dofs array

Args: - m: poloidal mode number - n: toroidal mode number (normalized by Nfp) - mpol: resolution of dofs array. If None (by default), use self.mpol - ntor: resolution of dofs array. If None (by default), use self.ntor - even: set to True to get vnc. Default is False

get_vnc(m, n)
get_vns(m, n)
set_vnc(m, n, value)
set_vns(m, n, value)

simsopt.field.sampling module

simsopt.field.sampling.draw_uniform_on_curve(curve, nsamples, safetyfactor=10)

Uses rejection sampling to sample points on a curve. Warning: assumes that the underlying quadrature points on the Curve are uniformly distributed.

Parameters:
  • curve – The simsopt.geo.curve.Curve to spawn the particles on.

  • nsamples – number of samples.

  • safetyfactor – how many more samples than nsamples to generate for rejection/acceptance.

simsopt.field.sampling.draw_uniform_on_surface(surface, nsamples, safetyfactor=10)

Uses rejection sampling to sample points on a surface. Warning: assumes that the underlying quadrature points on the surface are uniformly distributed.

Parameters:
  • surface – The simsopt.geo.surface.Surface to spawn the particles on.

  • nsamples – number of samples.

  • safetyfactor – how many more samples than nsamples to generate for rejection/acceptance.

simsopt.field.tracing module